Library
Public API
Creating an Assembly
GXBeam.curve_length
— Functioncurve_length(start, stop, curvature)
Calculate the length of a curve given its endpoints and its curvature vector.
GXBeam.discretize_beam
— Functiondiscretize_beam(L, start, discretization; frame, curvature)
Discretize a beam of length L
located at start
according to the discretization provided in discretization
Return the lengths, endpoints, midpoints, and reference frame of the beam elements.
Arguments
L
: Beam lengthstart
: Beam starting pointdiscretization
: Number of beam elements, or the normalized endpoints of each beam element, with values ranging from 0 to 1.
Keyword Arguments
frame
: Reference frame at the start of the beam element, represented by a 3x3 transformation matrix from the undeformed local frame to the body frame.curvature
: curvature vector
discretize_beam(start, stop, discretization; frame, curvature)
Discretize a beam from start
to stop
according to the discretization provided in discretization
.
Return the lengths, endpoints, midpoints, and reference frame of the beam elements.
Arguments
start
: Beam starting pointstop
: Beam ending pointdiscretization
: Number of beam elements, or the normalized endpoints of each beam element, with values ranging from 0 to 1.
Keyword Arguments
frame
: Reference frame at the start of the beam element, represented by a 3x3 transformation matrix from the undeformed local frame to the body frame.curvature
: curvature vector
GXBeam.Assembly
— TypeAssembly{TF, TP<:AbstractVector{<:AbstractVector{TF}},
TC<:AbstractVector{<:Integer}, TE<:AbstractVector{Element{TF}}}
Composite type that defines an assembly of connected nonlinear beam elements.
Fields
points
: Array of all beam element endpointsstart
: Array containing point index where each beam element startsstop
: Array containing point index where each beam element stopselements
: Array containing beam element definitions (seeElement
)
GXBeam.Assembly
— MethodAssembly(points, start, stop; kwargs...)
Construct an assembly of connected nonlinear beam elements. Beam lengths and midpoints may be manually specified in case beam elements are curved rather than straight.
Arguments
points
: Array of all beam element endpointsstart
: Array containing point indices where each beam element startsstop
: Array containing point indices where each beam element stops
Keyword Arguments
frames
: Array of (3 x 3) tranformation matrices for each beam element. Transforms from the local undeformed beam frame to the global frame. Defaults to the identity matrix.compliance
: Array of (6 x 6) compliance matrices for each beam element. Defaults tozeros(6,6)
for each beam element.mass
: Array of (6 x 6) mass matrices for each beam element. Defaults tozeros(6,6)
for each beam element.damping
: Array of (6) structural damping coefficients for each beam element. Defaults tofill(0.01, 6)
for each beam element.lengths
: Array containing the length of each beam element. Defaults to the distance between beam endpoints.midpoints
: Array containing the midpoint of each beam element. Defaults to the average of the beam element endpoints.
Section Properties
GXBeam.Material
— TypeMaterial(E1, E2, E3, G12, G13, G23, nu12, nu13, nu23, rho,
S1t=1.0, S1c=1.0, S2t=1.0, S2c=1.0, S3t=1.0, S3c=1.0, S12=1.0, S13=1.0, S23=1.0)
Orthotropic material properties. Axis 1 is the main ply axis, axis 2 is the transverse ply axis, and axis 3 is normal to the ply. For a fiber orientation of zero, axis 1 is along the beam axis.
Arguments
Ei::float
: Young's modulus along ith axis.Gij::float
: Shear modulinuij::float
: Poisson's ratio. $nu_ij E_j = nu_ji E_i$rho::float
: DensitySit::float
: Tensile strength in ith directionSic::float
: Compressive strength in ith directionSij::float
: Strength in ij direction
GXBeam.Node
— TypeNode(x, y)
A node in the finite element mesh at location x, y. If assembled in a vector, the vector index is the node number.
Arguments
x::float
: x-location of node in global coordinate systemy::float
: y-location of node in global coordinate system
GXBeam.MeshElement
— TypeMeshElement(nodenum, material, theta)
An element in the mesh, consisting of four ordered nodes, a material, and a fiber orientation.
Arguments
nodenum::Vector{integer}
: Node indices, defined counterclockwise starting from the bottom left node (as defined in the local coordinate system, see the figure).material::Material
: Material propertiestheta::float
: Fiber orientation
GXBeam.Layer
— TypeLayer(material, t, theta)
A layer (could be one ply or many plys of same material). A layup is a vector of layers.
Arguments
material::Material
: material of layert::float
: thickness of layertheta::float
: fiber orientation (rad)
GXBeam.afmesh
— Functionafmesh(xaf, yaf, chord, twist, paxis, xbreak, webloc, segments, webs; ds=nothing, dt=nothing, ns=nothing, nt=nothing, wns=4, wnt=nothing)
Create structural mesh for airfoil. The airfoil coordinates define the meshing density tangential to the airfoil. Whereas the number of layers defines the resolution normal to the airfoil. All segments are meshed with the same resolution in the normal direction, using the number of grid points as defined by segment with the most layers.
Arguments
xaf, yaf::Vector{float}
: points defining airfoil start at trailing edge and traverse counterclockwise back to trailing edge (can be blunt or sharp T.E.)chord::float
: chord lengthtwist::float
: twist angle (rad)paxis::float
: pitch axis (normalized by chord). e.g., 0.25 means twist about quarter chord.xbreak::Vector{float}
: x-locations, normalized by chord, defining break points between segments. must start at 0 and end at 1. e.g., [0, 0.2, 0.4, 0.7, 1.0] defines 4 segments.webloc::Vector{float}
: x-locations, normalized by chord, defining web centers (and length of vector is number of webs). e.g., [0.25, 0.55] means there is a web at 25% chord and a second web at 55% chord.segments::Vector{Vector{Layer}}
: A Layer defines a ply (or multiple plys with the same material/orientation). At a given x location the ply stack (segment) is defined a vector of layers starting from outside surface towards inside. Segments then is a vector of these segments that defines properties between segments as defined by xbreak.webs::Vector{Vector{Layer}}
: same structure as segments, except each inner vector is from left to right (although this is usually symmetric), and each outer vector is for a separate webds::float
: if provided, airfoil spacing will be resampling with approximately this spacing, normalized by chord. e.g., 0.01 will have points on the airfoil roughly 1% chord apart.dt::float
: if provided, thickness will be resampled with this maximum mesh size (thickness is absolute). Note that the total number of cells remains constant along airfoil, so most thicknesses will be much less. e.g., 0.01 will target a maximum mesh thickness of 0.01 (absolute).ns::vector{int}
: if provided, rather than use a targert size ds, we specify the number of cells to use in each segment. This is desirable for gradient-based optimization, if airfoil coordinates are changed, so that during resizing operations the mesh stretch/shrinks rather than experiencing discrete jumps. For example, ns=[15, 20, 40, 30] would use 15 elements between xbreak[1] and xbreak[2] and so on.nt::vector{vector{int}}
: if provided, defines how many elements to use across tangential direction. Again, prefered over dt for gradient-based optimization, if the thicknesses are changed during optimization. each entry defines how many cells to put in that layer following order of original layup. for example, nt=[[1, 2, 1], [1, 3]] would use 1 element, 2 elements (subdivide), then 1 elements over first sector, and so on.wns::int
: discretization level for number of elements vertically along web.wnt::vector{vector{int}}
: same definition as nt but for the webs
Returns
nodes::Vector{Node{Float64}}
: nodes for this meshelements::Vector{MeshElement{Float64}}
: elements for this mesh
GXBeam.initialize_cache
— Functioninitialize_cache([TF,] nodes, elements)
Construct the cache for the section property calculations.
GXBeam.compliance_matrix
— Functioncompliance_matrix(nodes, elements; cache=initialize_cache(nodes, elements),
gxbeam_order=true, shear_center=true)
Compute compliance matrix given a finite element mesh described by nodes and elements.
Arguments
nodes
: Vector containing all the nodes in the meshelements
: Vector containing all the elements in the meshcache::SectionCache
: A pre-allocated cache which may be passed in to reduce allocations across multiple calls to this function when the number of nodes, number of elements, and connectivity remain the same.gxbeam_order::Bool
: Indicates whether the compliance matrix should be provided in the order expected by GXBeam (rather than the internal ordering used by the section analysis)shear_center::Bool
: Indicates whether the compliance matrix should be provided about the shear center
Returns
S::Matrix
: compliance matrixsc::Vector{float}
: x, y location of shear center (location where a transverse/shear force will not produce any torsion, i.e., beam will not twist)tc::Vector{float}
: x, y location of tension center, aka elastic center, aka centroid (location where an axial force will not produce any bending, i.e., beam will remain straight)
GXBeam.mass_matrix
— Functionmass_matrix(nodes, elements)
Compute mass matrix for the structure using GXBeam ordering.
Returns
M::Matrix
: mass matrixmc::Vector{float}
: x, y location of mass center
GXBeam.plotmesh
— Functionplotmesh(nodes, elements, pyplot; plotnumbers=false)
plot nodes and elements for a quick visualization. Need to pass in a PyPlot object as PyPlot is not loaded by this package.
GXBeam.strain_recovery
— Functionstrain_recovery(F, M, nodes, elements, cache)
Compute stresses and strains at each element in cross section.
Arguments
F::Vector(3)
: force at this cross section in x, y, z directionsM::Vector(3)
: moment at this cross section in x, y, z directionsnodes::Vector{Node{TF}}
: all the nodes in the meshelements::Vector{MeshElement{TF}}
: all the elements in the meshcache::SectionCache
: needs to reuse data from the compliance solve (thus must initialize cache and pass it to both compliance and this function)
Returns
strain_b::Vector(6, ne)
: strains in beam coordinate system for each element. order: xx, yy, zz, xy, xz, yzstress_b::Vector(6, ne)
: stresses in beam coordinate system for each element. order: xx, yy, zz, xy, xz, yzstrain_p::Vector(6, ne)
: strains in ply coordinate system for each element. order: 11, 22, 33, 12, 13, 23stress_p::Vector(6, ne)
: stresses in ply coordinate system for each element. order: 11, 22, 33, 12, 13, 23
GXBeam.plotsoln
— Functionplotsoln(nodes, elements, soln, pyplot)
plot stress/strain on mesh soln could be any vector that is of length # of elements, e.g., sigma_b[3, :] Need to pass in a PyPlot object as PyPlot is not loaded by this package.
GXBeam.tsai_wu
— Functiontsai_wu(stress_p, elements)
Tsai Wu failure criteria
Arguments
stress_p::vector(6, ne)
: stresses in ply coordinate systemelements::Vector{MeshElement{TF}}
: all the elements in the mesh
Returns
failure::vector(ne)
: tsai-wu failure criteria for each element. fails if >= 1
Defining Point Masses
GXBeam.PointMass
— TypePointMass{T}
Type which contains the aggregated inertial properties of one or more point masses which are rigidly attached to the center of an element.
Fields:
mass
: Mass matrix corresponding to the point masses.
GXBeam.PointMass
— MethodPointMass(m, p, I)
Define a point mass given its mass m
, offset p
, and inertia matrix I
GXBeam.combine_masses
— Functioncombine_masses(masses)
Combine the point masses in the iterable collection masses
Defining Distributed Loads
GXBeam.DistributedLoads
— TypeDistributedLoads{T}
Type which contains pre-integrated distributed forces and moments applied to a beam element.
Fields
- f1: Integrated non-follower distributed force corresponding to the start of the beam element.
- f2: Integrated non-follower distributed force corresponding to the end of the beam element.
- m1: Integrated non-follower distributed moment corresponding to the start of the beam element.
- m2: Integrated non-follower distributed moment corresponding to the end of the beam element.
- f1_follower: Integrated follower distributed force corresponding to the start of the beam element.
- f2_follower: Integrated follower distributed force corresponding to the end of the beam element.
- m1_follower: Integrated follower distributed moment corresponding to the start of the beam element.
- m2_follower: Integrated follower distributed moment corresponding to the end of the beam element.
GXBeam.DistributedLoads
— MethodDistributedLoads(assembly, ielem; kwargs...)
Pre-integrate distributed loads on a beam element for use in an analysis.
Arguments
assembly
: Beam element assembly (of typeAssembly
)ielem
: Beam element index
Keyword Arguments
s1 = 0.0
: Start of the beam element (used solely for integrating the distributed loads)s2 = 1.0
: End of the beam element (used solely for integrating the distributed loads)method = (f, s1, s2) -> gauss_quadrature(f, s1, s2)
: Method which integrates function
f
from s1
to s2
. Defaults to the Gauss-Legendre quadrature with 4 points on each element.
fx = (s) -> 0.0
: Distributed x-direction forcefy = (s) -> 0.0
: Distributed y-direction forcefz = (s) -> 0.0
: Distributed z-direction forcemx = (s) -> 0.0
: Distributed x-direction momentmy = (s) -> 0.0
: Distributed y-direction momentmz = (s) -> 0.0
: Distributed z-direction momentfx_follower = (s) -> 0.0
: Distributed x-direction follower forcefy_follower = (s) -> 0.0
: Distributed y-direction follower forcefz_follower = (s) -> 0.0
: Distributed z-direction follower forcemx_follower = (s) -> 0.0
: Distributed x-direction follower momentmy_follower = (s) -> 0.0
: Distributed y-direction follower momentmz_follower = (s) -> 0.0
: Distributed z-direction follower moment
GXBeam.combine_loads
— Functioncombine_loads(loads)
Combine the distributed loads in the iterable collection loads
Defining Prescribed Conditions
GXBeam.PrescribedConditions
— TypePrescribedConditions{T}
Type which defines the prescribed displacements and loads at a point.
Fields:
pd
: Flag for each degree of freedom indicating whether displacements are prescribedpl
: Flag for each degree of freedom indicating whether loads are prescribedu
: Linear displacementtheta
: Angular displacementF
: External forcesM
: External momentsFf
: Follower forcesMf
: Follower moments
GXBeam.PrescribedConditions
— MethodPrescribedConditions(; kwargs...)
Define the prescribed conditions at a point. Individual prescribed conditions may be assigned as either a scalar parameter or as a function of time.
Prescribed Wiener-Milenkovic parameters must satisfy the following inequality: sqrt(thetax^2 + thetay^2 + theta_z^2) <= 4. Note that this restriction still allows all possible rotations to be represented.
Note that if displacements and loads corresponding to the same degree of freedom are prescribed at the same point, the global body-fixed acceleration corresponding to the same degree of freedom will be modified to attempt to satisfy both conditions.
Keyword Arguments
ux
: Prescribed x-displacement (in the body frame)uy
: Prescribed y-displacement (in the body frame)uz
: Prescribed z-displacement (in the body frame)theta_x
: Prescribed first Wiener-Milenkovic parametertheta_y
: Prescribed second Wiener-Milenkovic parametertheta_z
: Prescribed third Wiener-Milenkovic parameterFx
: Prescribed x-direction forceFy
: Prescribed y-direction forceFz
: Prescribed z-direction forceMx
: Prescribed x-direction momentMy
: Prescribed y-direction momentMz
: Prescribed z-direction momentFx_follower
: Prescribed x-direction follower forceFy_follower
: Prescribed y-direction follower forceFz_follower
: Prescribed z-direction follower forceMx_follower
: Prescribed x-direction follower momentMy_follower
: Prescribed y-direction follower momentMz_follower
: Prescribed z-direction follower moment
Pre-Initializing Memory for an Analysis
GXBeam.AbstractSystem
— TypeAbstractSystem
Supertype for types which contain the system state, residual vector, and jacobian matrix.
GXBeam.StaticSystem
— TypeStaticSystem{TF, TV<:AbstractVector{TF}, TM<:AbstractMatrix{TF}} <: AbstractSystem
Contains the system state, residual vector, and jacobian matrix for a static system.
GXBeam.StaticSystem
— MethodStaticSystem([TF=eltype(assembly),] assembly; kwargs...)
Initialize an object of type StaticSystem
.
Arguments:
TF:
(optional) Floating point type, defaults to the floating point type ofassembly
assembly
: Assembly of rigidly connected nonlinear beam elements
Keyword Arguments
force_scaling
: Factor used to scale system forces/moments internally. If not specified, a suitable default will be chosen based on the entries of the beam element compliance matrices.
GXBeam.DynamicSystem
— TypeDynamicSystem{TF, TV<:AbstractVector{TF}, TM<:AbstractMatrix{TF}} <: AbstractSystem
Contains the system state, residual vector, and jacobian matrix for a dynamic system.
GXBeam.DynamicSystem
— MethodDynamicSystem([TF=eltype(assembly),] assembly; kwargs...)
Initialize an object of type DynamicSystem
.
Arguments:
TF:
(optional) Floating point type, defaults to the floating point type ofassembly
assembly
: Assembly of rigidly connected nonlinear beam elements
Keyword Arguments
force_scaling
: Factor used to scale system forces/moments internally. If not specified, a suitable default will be chosen based on the entries of the beam element compliance matrices.
GXBeam.ExpandedSystem
— TypeExpandedSystem{TF, TV<:AbstractVector{TF}, TM<:AbstractMatrix{TF}} <: AbstractSystem
Contains the system state, residual vector, and jacobian matrix for a constant mass matrix system.
GXBeam.ExpandedSystem
— MethodExpandedSystem([TF=eltype(assembly),] assembly; kwargs...)
Initialize an object of type ExpandedSystem
.
Arguments:
TF:
(optional) Floating point type, defaults to the floating point type ofassembly
assembly
: Assembly of rigidly connected nonlinear beam elements
Keyword Arguments
force_scaling
: Factor used to scale system forces/moments internally. If not specified, a suitable default will be chosen based on the entries of the beam element compliance matrices.
GXBeam.reset_state!
— Functionreset_state!(system)
Sets the state variables in system
to zero.
GXBeam.set_state!
— Functionset_state!([x,] system::StaticSystem, prescribed_conditions; kwargs...)
Set the state variables in system
(or in the vector x
) to the provided values.
Keyword Arguments
u
: Vector containing the linear displacement of each point.theta
: Vector containing the angular displacement of each point.F
: Vector containing the externally applied forces acting on each pointM
: Vector containing the externally applied moments acting on each pointFi
: Vector containing internal forces for each beam elementMi
: Vector containing internal moments for each beam element
set_state!([x,] system::DynamicSystem, prescribed_conditions; kwargs...)
Set the state variables in system
(or in the vector x
) to the provided values.
Keyword Arguments
u
: Vector containing the linear displacement of each point.theta
: Vector containing the angular displacement of each point.V
: Vector containing the linear velocity of each point.Omega
Vector containing the angular velocity of each pointF
: Vector containing the externally applied forces acting on each pointM
: Vector containing the externally applied moments acting on each pointFi
: Vector containing internal forces for each beam element (in the deformed element frame)Mi
: Vector containing internal moments for each beam element (in the deformed element frame)
set_state!([x,] system::ExpandedSystem, prescribed_conditions; kwargs...)
Set the state variables in system
(or in the vector x
) to the provided values.
Keyword Arguments
u
: Vector containing the linear displacement of each point.theta
: Vector containing the angular displacement of each point.V
: Vector containing the linear velocity of each point in the deformed point frameOmega
Vector containing the angular velocity of each point in the deformed point frameF
: Vector containing the externally applied forces acting on each pointM
: Vector containing the externally applied moments acting on each pointF1
: Vector containing resultant forces at the start of each beam element (in the deformed element frame)M1
: Vector containing resultant moments at the start of each beam element (in the deformed element frame)F2
: Vector containing resultant forces at the end of each beam element (in the deformed element frame)M2
: Vector containing resultant moments at the end of each beam element (in the deformed element frame)V_e
: Vector containing the linear velocity of each beam element in the deformed beam element reference frame.Omega_e
Vector containing the angular velocity of each beam element in the deformed beam element reference frame.
set_state!([x=system.x,] system, assembly, state; kwargs...)
Set the state variables in x
to the values in state
GXBeam.set_rate!
— Functionset_rate!([x=system.dx,] system, assembly, state::AssemblyState; kwargs...)
Set the state variable rates in dx
to the values in state
GXBeam.set_linear_displacement!
— Functionset_linear_displacement!([x,] system, prescribed_conditions, u, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the linear deflection of point ipoint
to the provided values.
GXBeam.set_angular_displacement!
— Functionset_angular_displacement!([x,] system, prescribed_conditions, theta, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the angular deflection of point ipoint
to the provided values.
GXBeam.set_external_forces!
— Functionset_external_forces!([x,] system, prescribed_conditions, F, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the external forces applied at point ipoint
to the provided values.
GXBeam.set_external_moments!
— Functionset_external_moments!([x,] system, prescribed_conditions, M, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the external moments applied at point ipoint
to the provided values.
GXBeam.set_linear_velocity!
— Functionset_linear_velocity!([x,] system, V, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the linear velocity of point ipoint
to the provided values.
GXBeam.set_angular_velocity!
— Functionset_angular_velocity!([x,] system, Omega, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the angular velocity of point ipoint
to the provided values.
GXBeam.set_internal_forces!
— Functionset_internal_forces!([x,] system::Union{StaticSystem,DynamicSystem}, Fi, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the internal forces of element ielem
to the provided values.
GXBeam.set_internal_moments!
— Functionset_internal_moments!([x,] system::Union{StaticSystem, DynamicSystem}, Mi, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the internal moments of element ielem
to the provided values.
GXBeam.set_start_forces!
— Functionset_start_forces!([x,] system, F1, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the resultant forces at the start of element ielem
to the provided values.
GXBeam.set_start_moments!
— Functionset_start_moments!([x,] system, M1, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the resultant moments at the start of element ielem
to the provided values.
GXBeam.set_end_forces!
— Functionset_end_forces!([x,] system, F2, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the resultant forces at the end of element ielem
to the provided values.
GXBeam.set_end_moments!
— Functionset_end_moments!([x,] system, M2, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the resultant moments at the end of element ielem
to the provided values.
GXBeam.set_point_linear_velocity!
— Functionset_point_linear_velocity!([x,] system::ExpandedSystem, V, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the linear velocity of point ipoint
to the provided values.
GXBeam.set_point_angular_velocity!
— Functionset_point_angular_velocity!([x,] system::ExpandedSystem, Omega, ipoint)
Set the state variables in system
(or in the vector x
) corresponding to the angular velocity of point ipoint
to the provided values.
GXBeam.set_element_linear_velocity!
— Functionset_element_linear_velocity!([x,] system::ExpandedSystem, V, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the linear velocity of beam element ielem
to the provided values.
GXBeam.set_element_angular_velocity!
— Functionset_element_angular_velocity!([x,] system::ExpandedSystem, Omega, ielem)
Set the state variables in system
(or in the vector x
) corresponding to the angular velocity of element ielem
to the provided values.
Performing an Analysis
GXBeam.static_analysis
— Functionstatic_analysis(assembly; kwargs...)
Perform a static analysis for the system of nonlinear beams contained in assembly
. Return the resulting system, the post-processed solution, and a convergence flag indicating whether the iteration procedure converged.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.gravity = [0,0,0]
: Gravity vector. If time varying, this input may be provided as a function of time.time = 0.0
: Current time or vector of times corresponding to each step. May be used in conjunction with time varying prescribed conditions, distributed loads, and body frame motion to gradually increase displacements and loads.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
which contains the initial state variables. If not provided (or set tonothing
), then the state variables stored insystem
(which default to zeros) will be used as the initial state variables.linear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planeshow_trace = false
: Flag indicating whether to display the solution progress.
Linear Analysis Keyword Arguments
update_linearization = false
: Flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables. Iffalse
, then the initial set of state variables will be used for the linearization.
Nonlinear Analysis Keyword Arguments (see NLsolve.jl)
method = :newton
: Solution method for nonlinear systems of equationslinesearch = LineSearches.BackTracking(maxstep=1e6)
: Line search for solving nonlinear systems of equationsftol = 1e-9
: Tolerance for solving the nonlinear system of equationsiterations = 1000
: Iteration limit when solving nonlinear systems of equations
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
, andgravity
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
GXBeam.static_analysis!
— Functionstatic_analysis!(system, assembly; kwargs...)
Pre-allocated version of static_analysis
.
GXBeam.steady_state_analysis
— Functionsteady_state_analysis(assembly; kwargs...)
Perform a steady-state analysis for the system of nonlinear beams contained in assembly
. Return the resulting system and a flag indicating whether the iteration procedure converged.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.linear_velocity = zeros(3)
: Prescribed linear velocity of the body frame. If time varying, this input may be provided as a function of time.angular_velocity = zeros(3)
: Prescribed angular velocity of the body frame. If time varying, this input may be provided as a function of time.linear_acceleration = zeros(3)
: Prescribed linear acceleration of the body frame. If time varying, this input may be provided as a function of time.angular_acceleration = zeros(3)
: Prescribed angular acceleration of the body frame. If time varying, this input may be provided as a function of time.gravity = [0,0,0]
: Gravity vector in the body frame. If time varying, this input may be provided as a function of time.time = 0.0
: Current time or vector of times corresponding to each step. May be used in conjunction with time varying prescribed conditions, distributed loads, and body frame motion to gradually increase displacements and loads.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
which contains the initial state variables. If not provided (or set tonothing
), then the state variables stored insystem
(which default to zeros) will be used as the initial state variables.structural_damping = false
: Indicates whether to enable structural dampinglinear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planeshow_trace = false
: Flag indicating whether to display the solution progress.
Linear Analysis Keyword Arguments
update_linearization = false
: Flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables. Iffalse
, then the initial set of state variables will be used for the linearization.
Nonlinear Analysis Keyword Arguments
method = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve the nonlinear system of equationsftol = 1e-9
: tolerance for solving the nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equations=
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
,linear_velocity
,angular_velocity
,linear_acceleration
,angular_acceleration
, andgravity
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
GXBeam.steady_state_analysis!
— Functionsteady_state_analysis!(system, assembly; kwargs...)
Pre-allocated version of steady_state_analysis
.
GXBeam.eigenvalue_analysis
— Functioneigenvalue_analysis(assembly; kwargs...)
Compute the eigenvalues and eigenvectors of the system of nonlinear beams contained in assembly
. Return the modified system, eigenvalues, eigenvectors, and a convergence flag indicating whether the corresponding steady-state analysis converged.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.linear_velocity = zeros(3)
: Prescribed linear velocity of the body frame. If time varying, this input may be provided as a function of time.angular_velocity = zeros(3)
: Prescribed angular velocity of the body frame. If time varying, this input may be provided as a function of time.linear_acceleration = zeros(3)
: Prescribed linear acceleration of the body frame. If time varying, this input may be provided as a function of time.angular_acceleration = zeros(3)
: Prescribed angular acceleration of the body frame. If time varying, this input may be provided as a function of time.gravity = [0,0,0]
: Gravity vector in the inertial frame. If time varying, this input may be provided as a function of time.time = 0.0
: Current time or vector of times corresponding to each step. May be used in conjunction with time varying prescribed conditions, distributed loads, and body frame motion to gradually increase displacements and loads.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
which contains the initial state variables. If not provided (or set tonothing
), then the state variables stored insystem
(which default to zeros) will be used as the initial state variables.structural_damping = false
: Indicates whether to enable structural dampinglinear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planeshow_trace = false
: Flag indicating whether to display the solution progress.
Linear Analysis Keyword Arguments
update_linearization = false
: Flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables. Iffalse
, then the initial set of state variables will be used for the linearization.
Nonlinear Analysis Keyword Arguments
method = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve the nonlinear system of equationsftol = 1e-9
: tolerance for solving the nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equations
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
,linear_velocity
,angular_velocity
,linear_acceleration
,angular_acceleration
, andgravity
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
Eigenvalue Analysis Keyword Arguments
nev = 6
: Number of eigenvalues to computesteady = reset_state && !linear
: Flag indicating whether the steady state solution should be found prior to performing the eigenvalue analysis.left = false
: Flag indicating whether to return left and right eigenvectors rather than just right eigenvectors.Uprev = nothing
: Previous left eigenvector matrix. May be provided in order to reorder eigenvalues based on results from a previous iteration.
GXBeam.eigenvalue_analysis!
— Functioneigenvalue_analysis!(system, assembly; kwargs...)
Pre-allocated version of eigenvalue_analysis
.
GXBeam.initial_condition_analysis
— Functioninitial_condition_analysis(assembly, t0; kwargs...)
Perform an analysis to obtain a consistent set of initial conditions. Return the resulting system and a flag indicating whether the iteration procedure converged.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.linear_velocity = zeros(3)
: Initial linear velocity of the body frame.angular_velocity = zeros(3)
: Initial angular velocity of the body frame.linear_acceleration = zeros(3)
: Initial linear acceleration of the body frame.angular_acceleration = zeros(3)
: Initial angular acceleration of the body frame.gravity = [0,0,0]
: Gravity vector in the inertial frame.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
which contains the initial state variables. If not provided (or set tonothing
), then the state variables stored insystem
(which default to zeros) will be used as the initial state variables.structural_damping = true
: Indicates whether to enable structural dampinglinear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planesteady_state=false
: Flag indicating whether to initialize by performing a steady state analysis.show_trace = false
: Flag indicating whether to display the solution progress.
Initial Condition Analysis Keyword Arguments
u0 = fill(zeros(3), length(assembly.points))
: Initial linear displacement of each point relative to the body frametheta0 = fill(zeros(3), length(assembly.points))
: Initial angular displacement of each point relative to the body frame (using Wiener-Milenkovic Parameters)V0 = fill(zeros(3), length(assembly.points))
: Initial linear velocity of each point relative to the body frameOmega0 = fill(zeros(3), length(assembly.points))
: Initial angular velocity of each point relative to the body frameVdot0 = fill(zeros(3), length(assembly.points))
: Initial linear acceleration of each point relative to the body frameOmegadot0 = fill(zeros(3), length(assembly.points))
: Initial angular acceleration of each point relative to the body frame
Linear Analysis Keyword Arguments
update_linearization = false
: Flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables. Iffalse
, then the initial set of state variables will be used for the linearization.
Nonlinear Analysis Keyword Arguments
method = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve the nonlinear system of equationsftol = 1e-9
: tolerance for solving the nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equations
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
,linear_velocity
,angular_velocity
,linear_acceleration
,angular_acceleration
,gravity
,u0
,theta0
,V0
,Omega0
,Vdot0
, andOmegadot0
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
GXBeam.initial_condition_analysis!
— Functioninitial_condition_analysis!(system, assembly, t0; kwargs...)
Pre-allocated version of initial_condition_analysis
.
GXBeam.time_domain_analysis
— Functiontime_domain_analysis(assembly, tvec; kwargs...)
Perform a time-domain analysis for the system of nonlinear beams contained in assembly
using the time vector tvec
. Return the final system, a post-processed solution history, and a convergence flag indicating whether the iteration procedure converged for every time step.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.linear_velocity = zeros(3)
: Prescribed linear velocity of the body frame.angular_velocity = zeros(3)
: Prescribed angular velocity of the body frame.linear_acceleration = zeros(3)
: Initial linear acceleration of the body frame.angular_acceleration = zeros(3)
: Initial angular acceleration of the body frame.gravity = [0,0,0]
: Gravity vector in the body frame. If time varying, this input may be provided as a function of time.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
, which defines the initial states and state rates corresponding to the analysis. By default, this input is calculated using eithersteady_state_analysis
orinitial_condition_analysis
.steady_state = false
: Flag indicating whether to compute the state variables corresponding to the keyword argumentinitial_state
usingsteady_state_analysis
(rather thaninitial_condition_analysis
).structural_damping = true
: Flag indicating whether to enable structural dampinglinear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planeshow_trace = false
: Flag indicating whether to display the solution progress.save = eachindex(tvec)
: Steps at which to save the time history
Initial Condition Analysis Arguments
u0 = fill(zeros(3), length(assembly.points))
: Initial linear displacement of each point in the body frametheta0 = fill(zeros(3), length(assembly.points))
: Initial angular displacement of each point in the body frame (using Wiener-Milenkovic Parameters)V0 = fill(zeros(3), length(assembly.points))
: Initial linear velocity of each point in the body frame excluding contributions from body frame motionOmega0 = fill(zeros(3), length(assembly.points))
: Initial angular velocity of each point in the body frame excluding contributions from body frame motionVdot0 = fill(zeros(3), length(assembly.points))
: Initial linear acceleration of each point in the body frame excluding contributions from body frame motionOmegadot0 = fill(zeros(3), length(assembly.points))
: Initial angular acceleration of each point in the body frame excluding contributions from body frame motion
Linear Analysis Keyword Arguments
update_linearization = false
: Flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables. Iffalse
, then the initial set of state variables will be used for the linearization.
Nonlinear Analysis Keyword Arguments
method = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve nonlinear systems of equationsftol = 1e-9
: tolerance for solving the nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equations
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
,linear_velocity
,angular_velocity
,linear_acceleration
,angular_acceleration
,gravity
,u0
,theta0
,V0
,Omega0
,Vdot0
, andOmegadot0
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
GXBeam.time_domain_analysis!
— Functiontime_domain_analysis!(system, assembly, tvec; kwargs...)
Pre-allocated version of time_domain_analysis
.
Post-Processing
GXBeam.AssemblyState
— TypeAssemblyState{TF, TP<:AbstractVector{PointState{TF}}, TE<:AbstractVector{ElementState{TF}}}
Struct for storing state variables for the points and elements in an assembly.
Fields:
points::TP
: Array ofPointState
for each point in the assemblyelements::TE
: Array ofElementState
for each element in the assembly
GXBeam.PointState
— TypePointState
Holds the state variables for a point
Fields:
u
: Linear displacementudot
: Linear displacement ratetheta
: Angular displacement (Wiener-Milenkovic parameters)thetadot
: Angular displacement rateV
: Linear velocityVdot
: Linear velocity rateOmega
: Angular velocityOmegadot
: Angular velocity rateF
: External forcesM
: External moments
GXBeam.ElementState
— TypeElementState
Holds the state variables for an element
Fields:
u
: Linear displacementudot
: Linear displacement ratetheta
: Angular displacement (Wiener-Milenkovic parameters)thetadot
: Angular displacement rateV
: Linear velocityVdot
: Linear velocity rateOmega
: Angular velocityOmegadot
: Angular velocity rateFi
: Internal forcesMi
: Internal moments
GXBeam.AssemblyState
— MethodAssemblyState(system, assembly; prescribed_conditions = Dict())
AssemblyState(x, system, assembly; prescribed_conditions = Dict())
AssemblyState(dx, x, system, assembly; prescribed_conditions = Dict())
Post-process the system state given the state vector x
and rate vector dx
. Return an object of type AssemblyState
that defines the state of the assembly for the time step.
If prescribed_conditions
is not provided, all point state variables are assumed to be displacements/rotations, rather than their actual identities as used in the analysis.
GXBeam.extract_element_state
— Functionextract_element_state(system, assembly, ielem; prescribed_conditions = Dict())
extract_element_state(x, system, assembly, ielem; prescribed_conditions = Dict())
extract_element_state(dx, x, system, assembly, ielem; prescribed_conditions = Dict())
Return the state variables corresponding to element ielem
(see ElementState
) given the state vector x
and rate vector dx
.
GXBeam.extract_element_states
— Functionextract_element_states(system, assembly, x = system.x, dx = system.dx;
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}())
Return the state variables corresponding to each element (see ElementState
) given the solution vector x
.
GXBeam.extract_point_state
— Functionextract_point_state(system, assembly, ipoint; prescribed_conditions = Dict())
extract_point_state(x, system, assembly, ipoint; prescribed_conditions = Dict())
extract_point_state(dx, x, system, assembly, ipoint; prescribed_conditions = Dict())
Return the state variables corresponding to point ipoint
(see PointState
) given the state vector x
and rate vector dx
.
If prescribed_conditions
is not provided, all point state variables are assumed to be displacements/rotations, rather than their actual identities as used in the analysis.
GXBeam.extract_point_states
— Functionextract_point_states(system, assembly; prescribed_conditions = Dict())
extract_point_states(x, system, assembly; prescribed_conditions = Dict())
extract_point_states(dx, x, system, assembly; prescribed_conditions = Dict())
Return the state variables corresponding to each point (see PointState
) given the state vector x
and rate vector dx
.
If prescribed_conditions
is not provided, all point state variables are assumed to be displacements/rotations, rather than their actual identities as used in the analysis.
GXBeam.linearize!
— Functionlinearize!(system, assembly; kwargs...)
Return the state variables, jacobian matrix, and mass matrix of a linearized system using the current system state vector.
General Keyword Arguments
prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}()
: A dictionary with keys corresponding to the points at which prescribed conditions are applied and values of typePrescribedConditions
which describe the prescribed conditions at those points. If time varying, this input may be provided as a function of time.distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: A dictionary with keys corresponding to the elements to which distributed loads are applied and values of typeDistributedLoads
which describe the distributed loads on those elements. If time varying, this input may be provided as a function of time.point_masses = Dict{Int,PointMass{Float64}}()
: A dictionary with keys corresponding to the points to which point masses are attached and values of typePointMass
which contain the properties of the attached point masses. If time varying, this input may be provided as a function of time.linear_velocity = zeros(3)
: Prescribed linear velocity of the body frame. If time varying, this input may be provided as a function of time.angular_velocity = zeros(3)
: Prescribed angular velocity of the body frame. If time varying, this input may be provided as a function of time.linear_acceleration = zeros(3)
: Prescribed linear acceleration of the body frame. If time varying, this input may be provided as a function of time.angular_acceleration = zeros(3)
: Prescribed angular acceleration of the body frame. If time varying, this input may be provided as a function of time.gravity = [0,0,0]
: Gravity vector in the inertial frame. If time varying, this input may be provided as a function of time.time = 0.0
: Current time or vector of times corresponding to each step. May be used in conjunction with time varying prescribed conditions, distributed loads, and body frame motion to gradually increase displacements and loads.
Control Flag Keyword Arguments
reset_state = true
: Flag indicating whether the system state variables should be set to zero prior to performing this analysis.initial_state = nothing
: Object of typeAssemblyState
which contains the initial state variables. If not provided (or set tonothing
), then the state variables stored insystem
(which default to zeros) will be used as the initial state variables.structural_damping = false
: Indicates whether to enable structural dampinglinear = false
: Flag indicating whether a linear analysis should be performed.two_dimensional = false
: Flag indicating whether to constrain results to the x-y planeshow_trace = false
: Flag indicating whether to display the solution progress.
Sensitivity Analysis Keyword Arguments
xpfunc = (x, p, t) -> (;)
: Similar topfunc
, except that parameters can also be defined as a function of GXBeam's state variables. Using this function forces the system jacobian to be computed using automatic differentiation and switches the nonlinear solver to a Newton-Krylov solver (with linesearch).pfunc = (p, t) -> (;)
: Function which returns a named tuple with fields corresponding to updated versions of the argumentsassembly
,prescribed_conditions
,distributed_loads
,point_masses
,linear_velocity
,angular_velocity
,linear_acceleration
,angular_acceleration
, andgravity
. Only fields contained in the resulting named tuple will be overwritten.p
: Sensitivity parameters, as defined in conjunction with the keyword argumentpfunc
. While not necessary, usingpfunc
andp
to define the arguments to this function allows automatic differentiation sensitivities to be computed more efficiently
GXBeam.solve_eigensystem
— Functionsolve_eigensystem(x, K, M, nev)
Return the eigenvalues and eigenvectors of a linearized system.
GXBeam.left_eigenvectors
— Functionleft_eigenvectors(system, λ, V)
left_eigenvectors(K, M, λ, V)
Compute the left eigenvector matrix U
for the system
using inverse power iteration given the eigenvalues λ
and the corresponding right eigenvector matrix V
.
The complex conjugate of each left eigenvector is stored in each row of the matrix U
Left and right eigenvectors satisfy the following M-orthogonality condition:
- u'Mv = 1 if u and v correspond to the same eigenvalue
- u'Mv = 0 if u and v correspond to different eigenvalues
This means that UMV = I
This function assumes that system
has not been modified since the eigenvalues and right eigenvectors were computed.
GXBeam.correlate_eigenmodes
— Functioncorrelate_eigenmodes(C)
Return the permutation and the associated corruption index vector which associates eigenmodes from the current iteration with those of the previous iteration given the correlation matrix C
.
The correlation matrix can take one of the following forms (in order of preference):
C = U_p*M*V
C = U*M_p*V_p
C = V_p'*V
C = V'*V_p
where U
is a matrix of conjugated left eigenvectors, M
is the system mass matrix, V
is a matrix of right eigenvectors, and ()_p
indicates a variable from the previous iteration.
Note that the following two forms of the correlation matrix seem to be significantly inferior to their counterparts listed above: C = U*M*V_p
and C = U_p*M_p*V
. This is likely due to the way in which the left eigenvector matrix is calculated.
The corruption index is the largest magnitude in a given row of C
that was not chosen divided by the magnitude of the chosen eigenmode. It is most meaningful when using one of the forms of the correlation matrix that uses left eigenvectors since correct eigenmodes will have magnitudes close to 1 and incorrect eigenmodes will have magnitudes close to 0.
If the new mode number is already assigned, the next highest unassigned mode number is used. In this case a corruption index higher than 1 will be returned, otherwise the values of the corruption index will always be bounded by 0 and 1.
See "New Mode Tracking Methods in Aeroelastic Analysis" by Eldred, Vankayya, and Anderson.
GXBeam.wiener_milenkovic
— Functionwiener_milenkovic(c)
Construct a Wiener-Milenkovic transformation matrix, given the three Wiener-Milenkovic parameters in c
.
Note that the corresponding rotation matrix is the transpose of this transformation matrix.
GXBeam.rotate
— Functionrotate(xyz, r, theta)
Rotate the vectors in xyz
about point r
using the Wiener-Milenkovic parameters in theta
.
GXBeam.rotate!
— Functionrotate!(xyz, r, theta)
Pre-allocated version of rotate
GXBeam.translate
— Functiontranslate(xyz, u)
Translate the points in xyz
by the displacements in u
.
GXBeam.translate!
— Functiontranslate!(xyz, u)
Pre-allocated version of translate
GXBeam.deform_cross_section
— Functiondeform_cross_section(xyz, r, u, theta)
Rotate the points in xyz
(of shape (3, :)) about point r
using the Wiener-Milenkovic parameters in theta
, then translate the points by the displacements in u
.
GXBeam.deform_cross_section!
— Functiondeform_cross_section!(xyz, r, u, theta)
Pre-allocated version of deform_cross_section
GXBeam.write_vtk
— Functionwrite_vtk(name, assembly::Assembly; kwargs...)
write_vtk(name, assembly::Assembly, state::AssemblyState; kwargs...)
write_vtk(name, assembly::Assembly, history::Vector{<:AssemblyState}], dt;
kwargs...)
Write the deformed geometry (and associated data) to a VTK file for visualization using ParaView.
The state
argument may be omitted to write the original geometry to a VTK file without any associated data.
If the solution time history
is provided, the time step must also be provided
Keyword Arguments
sections = nothing
: Cross section geometry corresponding to each point, defined in a frame aligned with the body frame but centered around the corresponding point. Defined as an array with shape(3, ncross, np)
wherencross
is the number of points in each cross section andnp
is the number of points.scaling=1.0
: Parameter to scale the deflections (only valid if state is provided)metadata=Dict()
: Dictionary of metadata for the file(s)
write_vtk(name, assembly::Assembly, [state::AssemblyState, ]λ::Number,
eigenstate::AssemblyState; scaling=1.0, mode_scaling=1.0, cycles=1,
steps=100)
Write a series of files corresponding to the elastic motion of the assembly
about the deformed state encoded in state
defined by the eigenvalue λ
and the eigenvector encoded in eigenstate
over the time period specified by time
.
The steady-state deflections can be scaled with scaling
and the eigenmode deflections can be scaled using mode_scaling
.
The current time is encoded in the metadata tag "time"
Private API
Math
GXBeam.tilde
— Functiontilde(x)
Construct the cross product operator matrix
GXBeam.transform_properties
— Functiontransform_properties(K, T)
Applies the transformation T
to the stiffness or mass matrix K
GXBeam.rotation_parameter_scaling
— Functionrotation_parameter_scaling(θ)
Extract a scaling parameter which may be multiplied by the angular parameters to yield the Wiener-Milenkovic rotation parameters. Use of this scaling parameter allows deflections greater than 360 degrees.
GXBeam.get_C
— Functionget_C(θ)
Returns the transformation matrix C
given the three angular parameters in θ
.
GXBeam.get_C_θ
— Functionget_C_θ([C, ] θ)
Calculate the derivative of the Wiener-Milenkovic transformation matrix C
with respect to each of the rotation parameters in θ
.
GXBeam.get_Q
— Functionget_Q(θ)
Calculate the matrix Q as defined in the paper "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu given the rotational parameters in θ
.
GXBeam.get_Q_θ
— Functionget_Q_θ(θ)
get_Q_θ(Q, θ)
Calculate the derivative of the matrix Q
with respect to each of the rotation parameters in θ
.
GXBeam.get_Qinv
— Functionget_Qinv(θ)
Calculate the matrix inverse Qinv
as defined in the paper "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu given the rotational parameters in θ
.
GXBeam.get_Qinv_θ
— Functionget_Qinv_θ(θ)
Calculate the derivative of the matrix inverse Qinv
with respect to each of the rotation parameters in θ
.
GXBeam.get_ΔQ
— Functionget_ΔQ(θ, Δθ [, Q])
Calculate the matrix ΔQ
for structural damping calculations
GXBeam.get_ΔQ_θ
— Functionget_ΔQ_θ(θ, Δθ, [Q, Q_θ1, Q_θ2, Q_θ3])
Calculate the derivative of the matrix ΔQ
with respect to each of the rotation parameters in θ
.
GXBeam.mul3
— Functionmul3(A_1, A_2, A_3, b)
Return the product of a 3x3x3 tensor represented by A_1
, A_2
, and A_3
with the vector b
.
GXBeam.gauss_quadrature
— Functiongauss_quadrature(f, a, b)
Default gauss-quadrature function used for integrating distributed loads.
Body Frame
GXBeam.update_body_acceleration_indices!
— Functionupdate_body_acceleration_indices!(system, prescribed_conditions)
update_body_acceleration_indices!(indices, prescribed_conditions)
Updates the state variable indices corresponding to the body frame accelerations to correspond to the provided prescribed conditions.
GXBeam.body_accelerations
— Functionbody_accelerations(system, x=system.x; linear_acceleration=zeros(3), angular_acceleration=zeros(3))
Extract the linear and angular acceleration of the body frame from the state vector, if applicable. Otherwise return the provided linear and angular acceleration. This function is applicable only for steady state and initial condition analyses.
Points
GXBeam.point_loads
— Functionpoint_loads(x, ipoint, icol, force_scaling, prescribed_conditions)
Extract the loads F
and M
of point ipoint
from the state variable vector or prescribed conditions.
GXBeam.point_load_jacobians
— Functionpoint_load_jacobians(x, ipoint, icol, force_scaling, prescribed_conditions)
Calculate the load jacobians F_θ
, F_F
, M_θ
, and M_M
of point ipoint
.
GXBeam.point_displacement
— Functionpoint_displacement(x, ipoint, icol_point, prescribed_conditions)
Extract the displacements u
and θ
of point ipoint
from the state variable vector or prescribed conditions.
GXBeam.point_displacement_jacobians
— Functionpoint_displacement_jacobians(ipoint, prescribed_conditions)
Calculate the displacement jacobians u_u
and θ_θ
of point ipoint
.
GXBeam.point_displacement_rates
— Functionpoint_displacement_rates(dx, ipoint, icol, prescribed_conditions)
Extract the displacement rates udot
and θdot
of point ipoint
from the rate variable vector.
GXBeam.point_velocities
— Functionpoint_velocities(x, ipoint, icol_point)
Extract the velocities V
and Ω
of point ipoint
from the state variable vector
GXBeam.initial_point_displacement
— Functioninitial_point_displacement(x, ipoint, icol_point, prescribed_conditions,
rate_vars)
Extract the displacements u
and θ
of point ipoint
from the state variable vector or prescribed conditions for an initial condition analysis.
GXBeam.initial_point_velocity_rates
— Functioninitial_point_velocity_rates(x, ipoint, icol_point, prescribed_conditions,
Vdot0, Ωdot0, rate_vars)
Extract the velocity rates Vdot
and Ωdot
of point ipoint
from the state variable vector or provided initial conditions. Note that Vdot
and Ωdot
in this case do not include any contributions resulting from body frame motion.
GXBeam.initial_point_displacement_jacobian
— Functioninitial_point_displacement_jacobian(ipoint, icol_point, prescribed_conditions,
rate_vars)
Extract the displacement jacobians u_u
and θ_θ
of point ipoint
from the state variable vector or prescribed conditions for an initial condition analysis.
GXBeam.initial_point_velocity_rate_jacobian
— Functioninitial_point_velocity_rate_jacobian(ipoint, icol_point, prescribed_conditions,
rate_vars)
Return the velocity rate jacobians Vdot_Vdot
and Ωdot_Ωdot
of point ipoint
. Note that Vdot
and Ωdot
in this case do not include any contributions resulting from body frame motion.
GXBeam.static_point_properties
— Functionstatic_point_properties(x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity)
Calculate/extract the point properties needed to construct the residual for a static analysis
GXBeam.steady_point_properties
— Functionsteady_point_properties(x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))
Calculate/extract the point properties needed to construct the residual for a steady state analysis
GXBeam.initial_point_properties
— Functioninitial_point_properties(x, indices, rate_vars, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate/extract the point properties needed to construct the residual for a time domain analysis initialization.
GXBeam.newmark_point_properties
— Functionnewmark_point_properties(x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
udot_init, θdot_init, Vdot_init, Ωdot_init, dt)
Calculate/extract the point properties needed to construct the residual for a newmark-scheme time stepping analysis
GXBeam.dynamic_point_properties
— Functiondynamic_point_properties(dx, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)
Calculate/extract the point properties needed to construct the residual for a dynamic analysis
GXBeam.expanded_steady_point_properties
— Functionexpanded_steady_point_properties(x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))
Calculate/extract the point properties needed to construct the residual for a constant mass matrix system.
GXBeam.expanded_dynamic_point_properties
— Functionexpanded_dynamic_point_properties(dx, x, indices, force_scaling, assembly,
ipoint, prescribed_conditions, point_masses, gravity, linear_velocity,
angular_velocity)
Calculate/extract the point properties needed to construct the residual for a constant mass matrix system.
GXBeam.point_velocity_residuals
— Functionpoint_velocity_residuals(properties)
Calculate the velocity residuals rV
and rΩ
for a point for a steady state analysis.
GXBeam.expanded_point_velocity_residuals
— Functionexpanded_point_velocity_residuals(properties)
Calculate the velocity residuals rV
and rΩ
for a point for a constant mass matrix system.
GXBeam.static_point_resultants
— Functionstatic_point_resultants(properties)
Calculate the net loads F
and M
applied at a point for a static analysis.
GXBeam.dynamic_point_resultants
— Functiondynamic_point_resultants(properties)
Calculate the net loads F
and M
applied at a point for a steady analysis.
GXBeam.expanded_point_resultants
— Functionexpanded_point_resultants(properties)
Calculate the net loads F
and M
applied at a point for a constant mass matrix system.
GXBeam.static_point_residual!
— Functionstatic_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity)
Calculate and insert the residual entries corresponding to a point for a static analysis into the system residual vector.
GXBeam.steady_point_residual!
— Functionsteady_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
linear_acceleration, angular_acceleration)
Calculate and insert the residual entries corresponding to a point for a steady state analysis into the system residual vector.
GXBeam.initial_point_residual!
— Functioninitial_point_residual!(resid, x, indices, rate_vars,
force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate and insert the residual entries corresponding to a point for the initialization of a time domain analysis into the system residual vector.
GXBeam.newmark_point_residual!
— Functionnewmark_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
udot_init, θdot_init, Vdot_init, Ωdot_init, dt)
Calculate and insert the residual entries corresponding to a point for a newmark-scheme time marching analysis into the system residual vector.
GXBeam.dynamic_point_residual!
— Functiondynamic_point_residual!(resid, dx, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)
Calculate and insert the residual entries corresponding to a point for a dynamic analysis into the system residual vector.
GXBeam.expanded_steady_point_residual!
— Functionexpanded_steady_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
linear_acceleration, angular_acceleration)
Calculate and insert the residual entries corresponding to a point into the system residual vector for a constant mass matrix system.
GXBeam.expanded_dynamic_point_residual!
— Functionexpanded_dynamic_point_residual!(resid, dx, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)
Calculate and insert the residual entries corresponding to a point into the system residual vector for a constant mass matrix system.
GXBeam.static_point_jacobian_properties
— Functionstatic_point_jacobian_properties(properties, x, indices, force_scaling, assembly,
ipoint, prescribed_conditions, point_masses, gravity)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a static analysis
GXBeam.steady_point_jacobian_properties
— Functionsteady_point_jacobian_properties(properties, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity,
ub_ub, θb_θb, vb_vb, ωb_ωb, ab_ab, αb_αb)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a steady state analysis
GXBeam.initial_point_jacobian_properties
— Functioninitial_point_jacobian_properties(properties, x, indices, rate_vars,
force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a Newmark scheme time marching analysis
GXBeam.newmark_point_jacobian_properties
— Functionnewmark_point_jacobian_properties(properties, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity,
udot_init, θdot_init, Vdot_init, Ωdot_init, dt)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a Newmark scheme time marching analysis
GXBeam.dynamic_point_jacobian_properties
— Functiondynamic_point_jacobian_properties(properties, dx, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a dynamic analysis
GXBeam.expanded_steady_point_jacobian_properties
— Functionexpanded_steady_point_jacobian_properties(properties, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a constant mass matrix system
GXBeam.expanded_dynamic_point_jacobian_properties
— Functionexpanded_dynamic_point_jacobian_properties(properties, dx, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity)
Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a constant mass matrix system
GXBeam.mass_matrix_point_jacobian_properties
— Functionmass_matrix_point_jacobian_properties(x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses)
Calculate/extract the point properties needed to calculate the mass matrix jacobian entries corresponding to a point
GXBeam.expanded_mass_matrix_point_jacobian_properties
— Functionexpanded_mass_matrix_point_jacobian_properties(assembly, ipoint, prescribed_conditions, point_masses)
Calculate/extract the point properties needed to calculate the mass matrix jacobian entries corresponding to a point for a constant mass matrix system
GXBeam.static_point_resultant_jacobians
— Functionstatic_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for a static analysis.
GXBeam.steady_point_resultant_jacobians
— Functionsteady_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for a steady state analysis.
GXBeam.initial_point_resultant_jacobians
— Functioninitial_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for the initialization of a time domain analysis.
GXBeam.dynamic_point_resultant_jacobians
— Functiondynamic_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for a dynamic analysis.
GXBeam.expanded_steady_point_resultant_jacobians
— Functionexpanded_steady_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for a constant mass matrix system
GXBeam.expanded_dynamic_point_resultant_jacobians
— Functionexpanded_dynamic_point_resultant_jacobians(properties)
Calculate the jacobians for the net loads F
and M
applied at a point for a constant mass matrix system
GXBeam.mass_matrix_point_resultant_jacobians
— Functionmass_matrix_point_resultant_jacobians(properties)
Calculate the mass matrix jacobians for the net loads F
and M
applied at a point
GXBeam.initial_point_velocity_jacobians
— Functioninitial_point_velocity_jacobians(properties)
Calculate the jacobians of the velocity residuals rV
and rΩ
of a point for the initialization of a time domain analysis.
GXBeam.newmark_point_velocity_jacobians
— Functionnewmark_point_velocity_jacobians(properties)
Calculate the jacobians of the velocity residuals rV
and rΩ
of a point for a Newmark scheme time-marching analysis.
GXBeam.dynamic_point_velocity_jacobians
— Functiondynamic_point_velocity_jacobians(properties)
Calculate the jacobians of the velocity residuals rV
and rΩ
of a point for a dynamic state analysis.
GXBeam.expanded_point_velocity_jacobians
— Functionexpanded_point_velocity_jacobians(properties)
Calculate the jacobians of the velocity residuals rV
and rΩ
of a point for a constant mass matrix system
GXBeam.mass_matrix_point_velocity_jacobians
— Functionmass_matrix_point_velocity_jacobians(properties)
Calculate the mass matrix jacobians of the velocity residuals rV
and rΩ
of a point
GXBeam.insert_static_point_jacobians!
— Functioninsert_static_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
resultants)
Insert the jacobian entries corresponding to a point for a static analysis into the system jacobian matrix.
GXBeam.insert_initial_point_jacobians!
— Functioninsert_initial_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
resultants, velocities)
Insert the jacobian entries corresponding to a point for the initialization of a time domain analysis into the system jacobian matrix.
GXBeam.insert_dynamic_point_jacobians!
— Functioninsert_dynamic_point_jacobians!(jacob, indices, force_scaling, ipoint,
properties, resultants, velocities)
Insert the jacobian entries corresponding to a point for a steady state analysis into the system jacobian matrix.
GXBeam.insert_expanded_steady_point_jacobians!
— Functioninsert_expanded_steady_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
resultants, velocities)
Insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix for a constant mass matrix system.
GXBeam.insert_expanded_dynamic_point_jacobians!
— Functioninsert_expanded_dynamic_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
resultants, velocities)
Insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix for a constant mass matrix system.
GXBeam.insert_mass_matrix_point_jacobians!
— Functioninsert_mass_matrix_point_jacobians!(jacob, gamma, indices, two_dimensional, force_scaling, ipoint,
properties, resultants, velocities)
Insert the mass matrix jacobian entries corresponding to a point into the system jacobian matrix.
GXBeam.static_point_jacobian!
— Functionstatic_point_jacobian!(jacob, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity)
Calculate and insert the jacobian entries corresponding to a point for a static analysis into the system jacobian matrix.
GXBeam.steady_point_jacobian!
— Functionsteady_point_jacobian!(jacob, x, indices, force_scaling, assembly,
ipoint, prescribed_conditions, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Calculate and insert the jacobian entries corresponding to a point for a steady state analysis into the system jacobian matrix.
GXBeam.initial_point_jacobian!
— Functioninitial_point_jacobian!(jacob, x, indices, rate_vars,
force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate and insert the jacobian entries corresponding to a point for the initialization of a time domain analysis into the system jacobian matrix.
GXBeam.newmark_point_jacobian!
— Functionnewmark_point_jacobian!(jacob, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
udot_init, θdot_init, Vdot_init, Ωdot_init, dt))
Calculate and insert the jacobian entries corresponding to a point for a Newmark scheme time marching analysis into the system jacobian matrix.
GXBeam.dynamic_point_jacobian!
— Functiondynamic_point_jacobian!(jacob, dx, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)
Calculate and insert the jacobian entries corresponding to a point for a dynamic analysis into the system jacobian matrix.
GXBeam.expanded_steady_point_jacobian!
— Functionexpanded_steady_point_jacobian!(jacob, x, indices, force_scaling,
assembly, ipoint, prescribed_conditions, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Calculate and insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix.
GXBeam.expanded_dynamic_point_jacobian!
— Functionexpanded_dynamic_point_jacobian!(jacob, dx, x, indices, force_scaling, assembly, ipoint,
prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)
Calculate and insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix.
GXBeam.mass_matrix_point_jacobian!
— Functionmass_matrix_point_jacobian!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
ipoint, prescribed_conditions)
Calculate and insert the mass_matrix jacobian entries corresponding to a point into the system jacobian matrix.
GXBeam.expanded_mass_matrix_point_jacobian!
— Functionexpanded_mass_matrix_point_jacobian!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
ipoint, prescribed_conditions, point_masses)
Calculate and insert the mass_matrix jacobian entries corresponding to a point into the system jacobian matrix for a constant mass matrix system
Elements
GXBeam.Element
— TypeElement{TF}
Composite type that defines a beam element's properties
Fields
L
: Beam element lengthx
: Beam element locationcompliance
: Beam element compliance matrixmass
: Beam element mass matrixCab
: Transformation matrix from the undeformed beam element frame to the body framemu
: Beam element damping coefficients
GXBeam.element_loads
— Functionelement_loads(x, ielem, icol_elem, force_scaling)
Extract the internal loads (F
, M
) for a beam element from the state variable vector. These loads are expressed in the deformed element frame.
GXBeam.expanded_element_loads
— Functionexpanded_element_loads(x, ielem, icol_elem, force_scaling)
Extract the internal loads of a beam element at its endpoints (F1
, M1
, F2
, M2
) from the state variable vector for a constant mass matrix system.
GXBeam.expanded_element_velocities
— Functionexpanded_element_velocities(x, ielem, icol_elem)
Extract the velocities of a beam element from the state variable vector for a constant mass matrix system.
GXBeam.static_element_properties
— Functionstatic_element_properties(x, indices, force_scaling, assembly, ielem,
prescribed_conditions, gravity)
Calculate/extract the element properties needed to construct the residual for a static analysis
GXBeam.steady_element_properties
— Functionsteady_element_properties(x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))
Calculate/extract the element properties needed to construct the residual for a steady state analysis
GXBeam.initial_element_properties
— Functioninitial_element_properties(x, indices, rate_vars, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate/extract the element properties needed to construct the residual for a time-domain analysis initialization
GXBeam.newmark_element_properties
— Functionnewmark_element_properties(x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
Vdot_init, Ωdot_init, dt)
Calculate/extract the element properties needed to construct the residual for a newmark- scheme time stepping analysis
GXBeam.dynamic_element_properties
— Functiondynamic_element_properties(dx, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity)
Calculate/extract the element properties needed to construct the residual for a dynamic analysis
GXBeam.expanded_steady_element_properties
— Functionexpanded_steady_element_properties(x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))
Calculate/extract the element properties needed to construct the residual for a constant mass matrix system
GXBeam.expanded_dynamic_element_properties
— Functionexpanded_dynamic_element_properties(dx, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity)
Calculate/extract the element properties needed to construct the residual for a constant mass matrix system
Missing docstring for GXBeam.compatability_residuals
. Check Documenter's build log for details.
Missing docstring for GXBeam.expanded_element_steady_velocity_residuals
. Check Documenter's build log for details.
Missing docstring for GXBeam.expanded_element_dynamic_velocity_residuals
. Check Documenter's build log for details.
GXBeam.expanded_element_equilibrium_residuals
— Functionexpanded_element_equilibrium_residuals(properties, distributed_loads, ielem)
Calculate the equilibrium residuals of a beam element for a constant mass matrix system.
GXBeam.static_element_resultants
— Functionstatic_element_resultants(properties, distributed_loads, ielem)
Calculate the resultant loads applied at each end of a beam element for a static analysis.
GXBeam.dynamic_element_resultants
— Functiondynamic_element_resultants(properties, distributed_loads, ielem)
Calculate the resultant loads applied at each end of a beam element for a steady state analysis.
GXBeam.expanded_element_resultants
— Functionexpanded_element_resultants(properties)
Calculate the resultant loads applied at each end of a beam element for a constant mass matrix system.
GXBeam.insert_element_residuals!
— Functioninsert_element_residuals!(resid, indices, force_scaling, assembly, ielem,
compatibility, resultants)
Insert the residual entries corresponding to a beam element into the system residual vector.
GXBeam.insert_expanded_element_residuals!
— Functioninsert_expanded_element_residuals!(resid, indices, force_scaling, assembly, ielem,
compatibility, velocities, equilibrium, resultants)
Insert the residual entries corresponding to a beam element into the system residual vector for a constant mass matrix system.
GXBeam.static_element_residual!
— Functionstatic_element_residual!(resid, x, indices, force_scaling, assembly, ielem,
prescribed_conditions, distributed_loads, gravity)
Calculate and insert the residual entries corresponding to a beam element for a static analysis into the system residual vector.
GXBeam.steady_element_residual!
— Functionsteady_element_residual!(resid, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, distributed_loads, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Calculate and insert the residual entries corresponding to a beam element for a steady state analysis into the system residual vector.
GXBeam.initial_element_residual!
— Functioninitial_element_residual!(resid, x, indices, rate_vars,
force_scaling, structural_damping, assembly, ielem, prescribed_conditions,
distributed_loads, gravity, linear_velocity, angular_velocity,
linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate and insert the residual entries corresponding to a beam element for the initialization of a time domain simulation into the system residual vector.
GXBeam.newmark_element_residual!
— Functionnewmark_element_residual!(resid, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, distributed_loads, gravity,
linear_velocity, angular_velocity, Vdot_init, Ωdot_init, dt)
Calculate and insert the residual entries corresponding to a beam element for a newmark-scheme time marching analysis into the system residual vector.
GXBeam.dynamic_element_residual!
— Functiondynamic_element_residual!(resid, dx, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, distributed_loads, gravity,
linear_velocity, angular_velocity)
Calculate and insert the residual entries corresponding to a beam element for a dynamic analysis into the system residual vector.
GXBeam.expanded_steady_element_residual!
— Functionexpanded_steady_element_residual!(resid, x, indices, force_scaling, structural_damping,
assembly, ielem, prescribed_conditions, distributed_loads, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Calculate and insert the residual entries corresponding to a beam element for a constant mass matrix system into the system residual vector.
GXBeam.expanded_dynamic_element_residual!
— Functionexpanded_dynamic_element_residual!(resid, dx, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity)
Calculate and insert the residual entries corresponding to a beam element for a constant mass matrix system into the system residual vector.
GXBeam.static_element_jacobian_properties
— Functionstatic_element_jacobian_properties(properties, x, indices, force_scaling,
assembly, ielem, prescribed_conditions, gravity)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a beam element for a static analysis
GXBeam.steady_element_jacobian_properties
— Functionsteady_element_jacobian_properties(properties, x, indices,
force_scaling, structural_damping, assembly, ielem, prescribed_conditions,
gravity)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a steady state analysis
GXBeam.initial_element_jacobian_properties
— Functioninitial_element_jacobian_properties(properties, x, indices, rate_vars,
force_scaling, structural_damping, assembly, ielem, prescribed_conditions, gravity,
u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a Newmark scheme time marching analysis
GXBeam.newmark_element_jacobian_properties
— Functionnewmark_element_jacobian_properties(properties, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, gravity,
Vdot_init, Ωdot_init, dt)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a Newmark scheme time marching analysis
GXBeam.dynamic_element_jacobian_properties
— Functiondynamic_element_jacobian_properties(properties, dx, x, indices,
force_scaling, assembly, ielem, prescribed_conditions, gravity)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a dynamic analysis
GXBeam.expanded_steady_element_jacobian_properties
— Functionexpanded_element_jacobian_properties(properties, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, gravity)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a steady state analysis
GXBeam.expanded_dynamic_element_jacobian_properties
— Functionexpanded_dynamic_element_jacobian_properties(properties, dx, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, gravity)
Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a dynamic analysis
GXBeam.mass_matrix_element_jacobian_properties
— Functionmass_matrix_element_jacobian_properties(x, indices, force_scaling,
assembly, ielem, prescribed_conditions)
Calculate/extract the element properties needed to calculate the mass matrix jacobian entries corresponding to a beam element
GXBeam.expanded_mass_matrix_element_jacobian_properties
— Functionexpanded_mass_matrix_element_jacobian_properties(assembly, ielem, prescribed_conditions)
Calculate/extract the element properties needed to calculate the mass matrix jacobian entries corresponding to a beam element for a constant mass matrix system
Missing docstring for GXBeam.expanded_element_steady_velocity_jacobians
. Check Documenter's build log for details.
Missing docstring for GXBeam.expanded_element_dynamic_velocity_jacobians
. Check Documenter's build log for details.
GXBeam.expanded_steady_element_equilibrium_jacobians
— Functionexpanded_steady_element_equilibrium_jacobians(properties)
Calculate the jacobians of the element equilibrium residuals for a constant mass matrix system.
GXBeam.expanded_dynamic_element_equilibrium_jacobians
— Functionexpanded_dynamic_element_equilibrium_jacobians(properties)
Calculate the jacobians of the element equilibrium residuals for a constant mass matrix system.
GXBeam.expanded_mass_matrix_element_equilibrium_jacobians
— Functionexpanded_mass_matrix_element_equilibrium_jacobians(properties)
Calculate the mass matrix jacobians for the resultant loads applied at each end of a beam element for a constant mass matrix system
GXBeam.static_element_resultant_jacobians
— Functionstatic_element_resultant_jacobians(properties, distributed_loads, ielem)
Calculate the jacobians for the resultant loads applied at each end of a beam element for a static analysis.
GXBeam.steady_element_resultant_jacobians
— Functionsteady_element_resultant_jacobians(properties, distributed_loads, ielem)
Calculate the jacobians for the resultant loads applied at each end of a beam element for a steady state analysis.
GXBeam.initial_element_resultant_jacobians
— Functioninitial_element_resultant_jacobians(properties, distributed_loads, ielem)
Calculate the jacobians for the resultant loads applied at each end of V beam element for the initialization of a time domain analysis.
GXBeam.dynamic_element_resultant_jacobians
— Functiondynamic_element_resultant_jacobians(properties, distributed_loads, ielem)
Calculate the jacobians for the resultant loads applied at each end of a beam element for a dynamic analysis.
GXBeam.expanded_element_resultant_jacobians
— Functionexpanded_element_resultant_jacobians(properties)
Calculate the jacobians for the resultant loads applied at each end of a beam element for a constant mass matrix system.
GXBeam.mass_matrix_element_resultant_jacobians
— Functionmass_matrix_element_resultant_jacobians(properties)
Calculate the mass matrix jacobians for the resultant loads applied at each end of a beam element
Missing docstring for GXBeam.expanded_mass_matrix_element_velocity_jacobians
. Check Documenter's build log for details.
GXBeam.static_element_jacobian!
— Functionstatic_element_jacobian!(jacob, x, indices, force_scaling,
assembly, ielem, prescribed_conditions, distributed_loads, gravity)
Calculate and insert the jacobian entries corresponding to a beam element for a static analysis into the system jacobian matrix.
GXBeam.steady_element_jacobian!
— Functionsteady_element_jacobian!(jacob, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Calculate and insert the jacobian entries corresponding to a beam element for a steady state analysis into the system jacobian matrix.
GXBeam.initial_element_jacobian!
— Functioninitial_element_jacobian!(jacob, x, indices, rate_vars, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, ub, θb, vb, ωb, ab, αb, u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Calculate and insert the jacobian entries corresponding to a beam element for the initialization of a time domain analysis into the system jacobian matrix.
GXBeam.newmark_element_jacobian!
— Functionnewmark_element_jacobian!(jacob, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity, Vdot_init, Ωdot_init, dt)
Calculate and insert the jacobian entries corresponding to a beam element for a Newmark-scheme time marching analysis into the system jacobian matrix.
GXBeam.dynamic_element_jacobian!
— Functiondynamic_element_jacobian!(jacob, dx, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity)
Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.
GXBeam.expanded_steady_element_jacobian!
— Functionexpanded_steady_element_jacobian!(jacob, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity, linear_acceleration,
angular_acceleration)
Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.
GXBeam.expanded_dynamic_element_jacobian!
— Functionexpanded_dynamic_element_jacobian!(jacob, dx, x, indices, force_scaling,
structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
gravity, linear_velocity, angular_velocity)
Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.
GXBeam.mass_matrix_element_jacobian!
— Functionmass_matrix_element_jacobian!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
ielem, prescribed_conditions)
Calculate and insert the mass_matrix jacobian entries corresponding to a beam element into the system jacobian matrix.
GXBeam.expanded_mass_matrix_element_jacobian!
— Functionexpanded_mass_matrix_element_jacobian!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
ielem, prescribed_conditions)
Calculate and insert the mass_matrix jacobian entries corresponding to a beam element into the system jacobian matrix for a constant mass matrix system
System
GXBeam.SystemIndices
— TypeSystemIndices
Structure for holding indices for accessing the state variables and equations associated with each point and beam element in a system.
GXBeam.default_force_scaling
— Functiondefault_force_scaling(assembly)
Defines a suitable default force scaling factor based on the nonzero elements of the compliance matrices in assembly
.
GXBeam.curve_triad
— Functioncurve_triad(Cab, k, s)
curve_triad(Cab, kkt, ktilde, kn, s)
Return the transformation matrix at s
along the length of the beam given the curvature vector k
and the initial transformation matrix Cab
.
GXBeam.curve_coordinates
— Functioncurve_coordiantes(r, Cab, k, s)
curve_coordinates(r, Cab, kkt, ktilde, kn, s)
Return the coordinates at s
along the length of the beam given the starting point r
, initial transformation matrix Cab
, and curvature vector k
.
GXBeam.static_system_residual!
— Functionstatic_system_residual!(resid, x, indices, two_dimensional, force_scaling,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity)
Populate the system residual vector resid
for a static analysis
GXBeam.initial_system_residual!
— Functioninitial_system_residual!(resid, x, indices, rate_vars1, rate_vars2,
two_dimensional, force_scaling, structural_damping, assembly, prescribed_conditions,
distributed_loads, point_masses, gravity, linear_velocity, angular_velocity,
linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Populate the system residual vector resid
for the initialization of a time domain simulation.
GXBeam.steady_system_residual!
— Functionsteady_system_residual!(resid, x, indices, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, vb, ωb, ab, αb)
Populate the system residual vector resid
for a steady state analysis
GXBeam.newmark_system_residual!
— Functionnewmark_system_residual!(resid, x, indices, two_dimensional, force_scaling, structural_damping,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
linear_velocity, angular_velocity, udot_init, θdot_init, Vdot_init, Ωdot_init, dt)
Populate the system residual vector resid
for a Newmark scheme time marching analysis.
GXBeam.dynamic_system_residual!
— Functiondynamic_system_residual!(resid, dx, x, indices, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, linear_velocity, angular_velocity)
Populate the system residual vector resid
for a general dynamic analysis.
GXBeam.expanded_steady_system_residual!
— Functionexpanded_steady_system_residual!(resid, x, indices, two_dimensional, force_scaling, structural_damping,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)
Populate the system residual vector resid
for a constant mass matrix system.
GXBeam.expanded_dynamic_system_residual!
— Functionexpanded_dynamic_system_residual!(resid, dx, x, indices, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, linear_velocity, angular_velocity)
Populate the system residual vector resid
for a constant mass matrix system.
GXBeam.static_system_jacobian!
— Functionstatic_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity)
Populate the system jacobian matrix jacob
for a static analysis
GXBeam.steady_system_jacobian!
— Functionsteady_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, linear_velocity, angular_velocity, linear_acceleration,
angular_acceleration)
Populate the system jacobian matrix jacob
for a steady-state analysis
GXBeam.initial_system_jacobian!
— Functioninitial_system_jacobian!(jacob, x, indices, rate_vars1, rate_vars2, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, linear_velocity, angular_velocity, linear_acceleration,
angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)
Populate the system jacobian matrix jacob
for the initialization of a time domain simulation.
GXBeam.newmark_system_jacobian!
— Functionnewmark_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling, structural_damping,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
linear_velocity, angular_velocity, udot_init, θdot_init, Vdot_init, Ωdot_init, dt)
Populate the system jacobian matrix jacob
for a Newmark scheme time marching analysis.
GXBeam.dynamic_system_jacobian!
— Functiondynamic_system_jacobian!(jacob, dx, x, indices, two_dimensional, force_scaling,
structural_damping, assembly, prescribed_conditions, distributed_loads,
point_masses, gravity, linear_velocity, angular_velocity)
Populate the system jacobian matrix jacob
for a general dynamic analysis.
GXBeam.expanded_steady_system_jacobian!
— Functionexpanded_steady_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling, structural_damping,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
ub_p, θb_p, vb_p, ωb_p, ab_p, αb_p)
Populate the system jacobian matrix jacob
for a general dynamic analysis with a constant mass matrix system.
GXBeam.expanded_dynamic_system_jacobian!
— Functionexpanded_dynamic_system_jacobian!(jacob, dx, x, indices, two_dimensional, force_scaling, structural_damping,
assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
linear_velocity, angular_velocity)
Populate the system jacobian matrix jacob
for a general dynamic analysis with a constant mass matrix system.
GXBeam.system_mass_matrix!
— Functionsystem_mass_matrix!(jacob, x, indices, two_dimensional, force_scaling, assembly,
prescribed_conditions, point_masses)
Calculate the jacobian of the residual expressions with respect to the state rates.
system_mass_matrix!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
prescribed_conditions, point_masses)
Calculate the jacobian of the residual expressions with respect to the state rates and add the result multiplied by gamma
to jacob
.
GXBeam.expanded_system_mass_matrix
— Functionexpanded_system_mass_matrix(system, assembly;
two_dimensional = false,
prescribed_conditions=Dict{Int, PrescribedConditions}(),
point_masses=Dict{Int, PointMass}())
Calculate the jacobian of the residual expressions with respect to the state rates for a constant mass matrix system.
GXBeam.expanded_system_mass_matrix!
— Functionexpanded_system_mass_matrix!(jacob, indices, two_dimensional, force_scaling, assembly, prescribed_conditions,
point_masses)
Calculate the jacobian of the residual expressions with respect to the state rates.
expanded_system_mass_matrix!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
prescribed_conditions, point_masses)
Calculate the jacobian of the residual expressions with respect to the state rates and add the result multiplied by gamma
to jacob
.
Index
GXBeam.AbstractSystem
GXBeam.Assembly
GXBeam.Assembly
GXBeam.AssemblyState
GXBeam.AssemblyState
GXBeam.DistributedLoads
GXBeam.DistributedLoads
GXBeam.DynamicSystem
GXBeam.DynamicSystem
GXBeam.Element
GXBeam.ElementState
GXBeam.ExpandedSystem
GXBeam.ExpandedSystem
GXBeam.Layer
GXBeam.Material
GXBeam.MeshElement
GXBeam.Node
GXBeam.PointMass
GXBeam.PointMass
GXBeam.PointState
GXBeam.PrescribedConditions
GXBeam.PrescribedConditions
GXBeam.StaticSystem
GXBeam.StaticSystem
GXBeam.SystemIndices
SciMLBase.DAEFunction
SciMLBase.DAEProblem
SciMLBase.ODEFunction
SciMLBase.ODEProblem
GXBeam.afmesh
GXBeam.body_accelerations
GXBeam.combine_loads
GXBeam.combine_masses
GXBeam.compliance_matrix
GXBeam.correlate_eigenmodes
GXBeam.curve_coordinates
GXBeam.curve_length
GXBeam.curve_triad
GXBeam.default_force_scaling
GXBeam.deform_cross_section
GXBeam.deform_cross_section!
GXBeam.discretize_beam
GXBeam.dynamic_element_jacobian!
GXBeam.dynamic_element_jacobian_properties
GXBeam.dynamic_element_properties
GXBeam.dynamic_element_residual!
GXBeam.dynamic_element_resultant_jacobians
GXBeam.dynamic_element_resultants
GXBeam.dynamic_point_jacobian!
GXBeam.dynamic_point_jacobian_properties
GXBeam.dynamic_point_properties
GXBeam.dynamic_point_residual!
GXBeam.dynamic_point_resultant_jacobians
GXBeam.dynamic_point_resultants
GXBeam.dynamic_point_velocity_jacobians
GXBeam.dynamic_system_jacobian!
GXBeam.dynamic_system_residual!
GXBeam.eigenvalue_analysis
GXBeam.eigenvalue_analysis!
GXBeam.element_loads
GXBeam.expanded_dynamic_element_equilibrium_jacobians
GXBeam.expanded_dynamic_element_jacobian!
GXBeam.expanded_dynamic_element_jacobian_properties
GXBeam.expanded_dynamic_element_properties
GXBeam.expanded_dynamic_element_residual!
GXBeam.expanded_dynamic_point_jacobian!
GXBeam.expanded_dynamic_point_jacobian_properties
GXBeam.expanded_dynamic_point_properties
GXBeam.expanded_dynamic_point_residual!
GXBeam.expanded_dynamic_point_resultant_jacobians
GXBeam.expanded_dynamic_system_jacobian!
GXBeam.expanded_dynamic_system_residual!
GXBeam.expanded_element_equilibrium_residuals
GXBeam.expanded_element_loads
GXBeam.expanded_element_resultant_jacobians
GXBeam.expanded_element_resultants
GXBeam.expanded_element_velocities
GXBeam.expanded_mass_matrix_element_equilibrium_jacobians
GXBeam.expanded_mass_matrix_element_jacobian!
GXBeam.expanded_mass_matrix_element_jacobian_properties
GXBeam.expanded_mass_matrix_point_jacobian!
GXBeam.expanded_mass_matrix_point_jacobian_properties
GXBeam.expanded_point_resultants
GXBeam.expanded_point_velocity_jacobians
GXBeam.expanded_point_velocity_residuals
GXBeam.expanded_steady_element_equilibrium_jacobians
GXBeam.expanded_steady_element_jacobian!
GXBeam.expanded_steady_element_jacobian_properties
GXBeam.expanded_steady_element_properties
GXBeam.expanded_steady_element_residual!
GXBeam.expanded_steady_point_jacobian!
GXBeam.expanded_steady_point_jacobian_properties
GXBeam.expanded_steady_point_properties
GXBeam.expanded_steady_point_residual!
GXBeam.expanded_steady_point_resultant_jacobians
GXBeam.expanded_steady_system_jacobian!
GXBeam.expanded_steady_system_residual!
GXBeam.expanded_system_mass_matrix
GXBeam.expanded_system_mass_matrix!
GXBeam.extract_element_state
GXBeam.extract_element_states
GXBeam.extract_point_state
GXBeam.extract_point_states
GXBeam.gauss_quadrature
GXBeam.get_C
GXBeam.get_C_θ
GXBeam.get_Q
GXBeam.get_Q_θ
GXBeam.get_Qinv
GXBeam.get_Qinv_θ
GXBeam.get_ΔQ
GXBeam.get_ΔQ_θ
GXBeam.initial_condition_analysis
GXBeam.initial_condition_analysis!
GXBeam.initial_element_jacobian!
GXBeam.initial_element_jacobian_properties
GXBeam.initial_element_properties
GXBeam.initial_element_residual!
GXBeam.initial_element_resultant_jacobians
GXBeam.initial_point_displacement
GXBeam.initial_point_displacement_jacobian
GXBeam.initial_point_jacobian!
GXBeam.initial_point_jacobian_properties
GXBeam.initial_point_properties
GXBeam.initial_point_residual!
GXBeam.initial_point_resultant_jacobians
GXBeam.initial_point_velocity_jacobians
GXBeam.initial_point_velocity_rate_jacobian
GXBeam.initial_point_velocity_rates
GXBeam.initial_system_jacobian!
GXBeam.initial_system_residual!
GXBeam.initialize_cache
GXBeam.insert_dynamic_point_jacobians!
GXBeam.insert_element_residuals!
GXBeam.insert_expanded_dynamic_point_jacobians!
GXBeam.insert_expanded_element_residuals!
GXBeam.insert_expanded_steady_point_jacobians!
GXBeam.insert_initial_point_jacobians!
GXBeam.insert_mass_matrix_point_jacobians!
GXBeam.insert_static_point_jacobians!
GXBeam.left_eigenvectors
GXBeam.linearize!
GXBeam.mass_matrix
GXBeam.mass_matrix_element_jacobian!
GXBeam.mass_matrix_element_jacobian_properties
GXBeam.mass_matrix_element_resultant_jacobians
GXBeam.mass_matrix_point_jacobian!
GXBeam.mass_matrix_point_jacobian_properties
GXBeam.mass_matrix_point_resultant_jacobians
GXBeam.mass_matrix_point_velocity_jacobians
GXBeam.mul3
GXBeam.newmark_element_jacobian!
GXBeam.newmark_element_jacobian_properties
GXBeam.newmark_element_properties
GXBeam.newmark_element_residual!
GXBeam.newmark_point_jacobian!
GXBeam.newmark_point_jacobian_properties
GXBeam.newmark_point_properties
GXBeam.newmark_point_residual!
GXBeam.newmark_point_velocity_jacobians
GXBeam.newmark_system_jacobian!
GXBeam.newmark_system_residual!
GXBeam.plotmesh
GXBeam.plotsoln
GXBeam.point_displacement
GXBeam.point_displacement_jacobians
GXBeam.point_displacement_rates
GXBeam.point_load_jacobians
GXBeam.point_loads
GXBeam.point_velocities
GXBeam.point_velocity_residuals
GXBeam.reset_state!
GXBeam.rotate
GXBeam.rotate!
GXBeam.rotation_parameter_scaling
GXBeam.set_angular_displacement!
GXBeam.set_angular_velocity!
GXBeam.set_element_angular_velocity!
GXBeam.set_element_linear_velocity!
GXBeam.set_end_forces!
GXBeam.set_end_moments!
GXBeam.set_external_forces!
GXBeam.set_external_moments!
GXBeam.set_internal_forces!
GXBeam.set_internal_moments!
GXBeam.set_linear_displacement!
GXBeam.set_linear_velocity!
GXBeam.set_point_angular_velocity!
GXBeam.set_point_linear_velocity!
GXBeam.set_rate!
GXBeam.set_start_forces!
GXBeam.set_start_moments!
GXBeam.set_state!
GXBeam.solve_eigensystem
GXBeam.static_analysis
GXBeam.static_analysis!
GXBeam.static_element_jacobian!
GXBeam.static_element_jacobian_properties
GXBeam.static_element_properties
GXBeam.static_element_residual!
GXBeam.static_element_resultant_jacobians
GXBeam.static_element_resultants
GXBeam.static_point_jacobian!
GXBeam.static_point_jacobian_properties
GXBeam.static_point_properties
GXBeam.static_point_residual!
GXBeam.static_point_resultant_jacobians
GXBeam.static_point_resultants
GXBeam.static_system_jacobian!
GXBeam.static_system_residual!
GXBeam.steady_element_jacobian!
GXBeam.steady_element_jacobian_properties
GXBeam.steady_element_properties
GXBeam.steady_element_residual!
GXBeam.steady_element_resultant_jacobians
GXBeam.steady_point_jacobian!
GXBeam.steady_point_jacobian_properties
GXBeam.steady_point_properties
GXBeam.steady_point_residual!
GXBeam.steady_point_resultant_jacobians
GXBeam.steady_state_analysis
GXBeam.steady_state_analysis!
GXBeam.steady_system_jacobian!
GXBeam.steady_system_residual!
GXBeam.strain_recovery
GXBeam.system_mass_matrix!
GXBeam.tilde
GXBeam.time_domain_analysis
GXBeam.time_domain_analysis!
GXBeam.transform_properties
GXBeam.translate
GXBeam.translate!
GXBeam.tsai_wu
GXBeam.update_body_acceleration_indices!
GXBeam.wiener_milenkovic
GXBeam.write_vtk